home *** CD-ROM | disk | FTP | other *** search
/ SGI Freeware 2002 November / SGI Freeware 2002 November - Disc 2.iso / dist / fw_gsl.idb / usr / freeware / info / gsl-ref.info-6.z / gsl-ref.info-6
Text File  |  2000-10-09  |  48KB  |  1,204 lines

  1. This is gsl-ref.info, produced by Makeinfo version 3.12h from
  2. gsl-ref.texi.
  3.  
  4. INFO-DIR-SECTION Scientific software
  5. START-INFO-DIR-ENTRY
  6. * gsl-ref: (gsl-ref).                   GNU Scientific Library - Reference
  7. END-INFO-DIR-ENTRY
  8.  
  9.    This file documents the GNU Scientific Library.
  10.  
  11.    Copyright (C) 1996, 1997, 1998, 1999 The GSL Project.
  12.  
  13.    Permission is granted to make and distribute verbatim copies of this
  14. manual provided the copyright notice and this permission notice are
  15. preserved on all copies.
  16.  
  17.    Permission is granted to copy and distribute modified versions of
  18. this manual under the conditions for verbatim copying, provided that
  19. the entire resulting derived work is distributed under the terms of a
  20. permission notice identical to this one.
  21.  
  22.    Permission is granted to copy and distribute translations of this
  23. manual into another language, under the above conditions for modified
  24. versions, except that this permission notice may be stated in a
  25. translation approved by the Foundation.
  26.  
  27. 
  28. File: gsl-ref.info,  Node: Exponential Integrals,  Next: Fermi-Dirac Function,  Prev: Exponential Function,  Up: Special Functions
  29.  
  30. Exponential Integrals
  31. =====================
  32.  
  33. Exponential Integrals
  34. ---------------------
  35.  
  36.    We have  E_1(x) := Re \int_1^\infty dt \exp(-xt)/t ;  E_2(x) := Re
  37. \int_1^\infty dt \exp(-xt)/t^2 .
  38.  
  39.  - Function: int gsl_sf_expint_E1_impl (double x, gsl_sf_result *
  40.           result)
  41.  - Function: int gsl_sf_expint_E1_e (double x, gsl_sf_result * result)
  42.      Domain: x != 0.0 Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW,
  43.      GSL_EUNDRFLW
  44.  
  45.  - Function: int gsl_sf_expint_E2_impl (double x, gsl_sf_result *
  46.           result)
  47.  - Function: int gsl_sf_expint_E2_e (double x, gsl_sf_result * result)
  48.      Domain: x != 0.0 Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW,
  49.      GSL_EUNDRFLW
  50.  
  51. Ei(x)
  52. -----
  53.  
  54.    We have  Ei(x) := PV \int_(-x)^\infty dt \exp(-t)/t .
  55.  
  56.  - Function: int gsl_sf_expint_Ei_impl (double x, gsl_sf_result *
  57.           result)
  58.  - Function: int gsl_sf_expint_Ei_e (double x, gsl_sf_result * result)
  59.      Domain: x != 0.0 Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW,
  60.      GSL_EUNDRFLW
  61.  
  62. Hyperbolic Integrals
  63. --------------------
  64.  
  65.    We have  Shi(x) = \int_0^x dt \sinh(t)/t ;  Chi(x) := Re[ M_EULER +
  66. log(x) + \int_0^x dt (Cosh[t]-1)/t .
  67.  
  68.  - Function: int gsl_sf_Shi_impl (double x, gsl_sf_result * result)
  69.  - Function: int gsl_sf_Shi_e (double x, gsl_sf_result * result)
  70.      Exceptional Return Values: GSL_EOVRFLW, GSL_EUNDRFLW
  71.  
  72.  - Function: int gsl_sf_Chi_impl (double x, gsl_sf_result * result)
  73.  - Function: int gsl_sf_Chi_e (double x, gsl_sf_result * result)
  74.      Domain: x != 0.0 Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW,
  75.      GSL_EUNDRFLW
  76.  
  77. Ei_3(x)
  78. -------
  79.  
  80.    We have  Ei_3(x) = \int_0^x dt \exp(-t^3) .
  81.  
  82.  - Function: int gsl_sf_expint_3_impl (double x, gsl_sf_result * result)
  83.  - Function: int gsl_sf_expint_3_e (double x, gsl_sf_result * result)
  84.      Domain: x >= 0.0 Exceptional Return Values: GSL_EDOM
  85.  
  86. Trigonometric Integrals
  87. -----------------------
  88.  
  89.    We have  Si(x) = \int_0^x dt \sin(t)/t ;  Ci(x) = -\int_x^\infty
  90. \cos(t)/t .
  91.  
  92.  - Function: int gsl_sf_Si_impl (const double x, gsl_sf_result * result)
  93.  - Function: int gsl_sf_Si_e (double x, gsl_sf_result * result)
  94.      Exceptional Return Values: none
  95.  
  96.  - Function: int gsl_sf_Ci_impl (const double x, gsl_sf_result * result)
  97.  - Function: int gsl_sf_Ci_e (double x, gsl_sf_result * result)
  98.      Domain: x > 0.0 Exceptional Return Values: GSL_EDOM
  99.  
  100. Arctangent Integral
  101. -------------------
  102.  
  103.    We have  AtanInt(x) = \int_0^x dt \arctan(t)/t .
  104.  
  105.  - Function: int gsl_sf_atanint_impl (double x, gsl_sf_result * result);
  106.  - Function: int gsl_sf_atanint_e (double x, gsl_sf_result * result);
  107.      Domain: Exceptional Return Values:
  108.  
  109. 
  110. File: gsl-ref.info,  Node: Fermi-Dirac Function,  Next: Gamma Function,  Prev: Exponential Integrals,  Up: Special Functions
  111.  
  112. Fermi-Dirac Function
  113. ====================
  114.  
  115. Complete Fermi-Dirac Integrals
  116. ------------------------------
  117.  
  118.    F_j(x)   := 1/Gamma[j+1] Integral[ t^j /(Exp[t-x] + 1),
  119. (t,0,Infinity)]
  120.  
  121.    /* Complete integral F_(-1)(x) = e^x / (1 + e^x)  *
  122.  
  123.  - Function: int gsl_sf_fermi_dirac_m1_impl (double x, gsl_sf_result *
  124.           result)
  125.  - Function: int gsl_sf_fermi_dirac_m1_e (double x, gsl_sf_result *
  126.           result)
  127.      Exceptional Return Values: GSL_EUNDRFLW
  128.  
  129.    /* Complete integral F_0(x) = ln(1 + e^x)  *
  130.  
  131.  - Function: int gsl_sf_fermi_dirac_0_impl (double x, gsl_sf_result *
  132.           result)
  133.  - Function: int gsl_sf_fermi_dirac_0_e (double x, gsl_sf_result *
  134.           result)
  135.      Exceptional Return Values: GSL_EUNDRFLW
  136.  
  137.    /* Complete integral F_1(x)  *
  138.  
  139.  - Function: int gsl_sf_fermi_dirac_1_impl (double x, gsl_sf_result *
  140.           result)
  141.  - Function: int gsl_sf_fermi_dirac_1_e (double x, gsl_sf_result *
  142.           result)
  143.      Exceptional Return Values: GSL_EUNDRFLW, GSL_EOVRFLW
  144.  
  145.    /* Complete integral F_2(x)  *
  146.  
  147.  - Function: int gsl_sf_fermi_dirac_2_impl (double x, gsl_sf_result *
  148.           result)
  149.  - Function: int gsl_sf_fermi_dirac_2_e (double x, gsl_sf_result *
  150.           result)
  151.      Exceptional Return Values: GSL_EUNDRFLW, GSL_EOVRFLW
  152.  
  153.    /* Complete integral F_j(x)  * for integer j  *
  154.  
  155.  - Function: int gsl_sf_fermi_dirac_int_impl (int j, double x,
  156.           gsl_sf_result * result)
  157.  - Function: int gsl_sf_fermi_dirac_int_e (int j, double x,
  158.           gsl_sf_result * result)
  159.      Exceptional Return Values: GSL_EUNDRFLW, GSL_EOVRFLW
  160.  
  161.    /* Complete integral F_(-1/2)(x)  *
  162.  
  163.  - Function: int gsl_sf_fermi_dirac_mhalf_impl (double x, gsl_sf_result
  164.           * result)
  165.  - Function: int gsl_sf_fermi_dirac_mhalf_e (double x, gsl_sf_result *
  166.           result)
  167.      Exceptional Return Values: GSL_EUNDRFLW, GSL_EOVRFLW
  168.  
  169.    /* Complete integral F_(1/2)(x)  *
  170.  
  171.  - Function: int gsl_sf_fermi_dirac_half_impl (double x, gsl_sf_result
  172.           * result)
  173.  - Function: int gsl_sf_fermi_dirac_half_e (double x, gsl_sf_result *
  174.           result)
  175.      Exceptional Return Values: GSL_EUNDRFLW, GSL_EOVRFLW
  176.  
  177.    /* Complete integral F_(3/2)(x)  *
  178.  
  179.  - Function: int gsl_sf_fermi_dirac_3half_impl (double x, gsl_sf_result
  180.           * result)
  181.  - Function: int gsl_sf_fermi_dirac_3half_e (double x, gsl_sf_result *
  182.           result)
  183.      Exceptional Return Values: GSL_EUNDRFLW, GSL_EOVRFLW
  184.  
  185. Incomplete Fermi-Dirac Integrals:
  186. ---------------------------------
  187.  
  188.    F_j(x,b) := 1/Gamma[j+1] Integral[ t^j /(Exp[t-x] + 1),
  189. (t,b,Infinity)]
  190.  
  191.    /* Incomplete integral F_0(x,b) = ln(1 + e^(b-x)) - (b-x)  *
  192.  
  193.  - Function: int gsl_sf_fermi_dirac_inc_0_impl (double x, double b,
  194.           gsl_sf_result * result)
  195.  - Function: int gsl_sf_fermi_dirac_inc_0_e (double x, double b,
  196.           gsl_sf_result * result)
  197.      Exceptional Return Values: GSL_EUNDRFLW, GSL_EDOM
  198.  
  199. 
  200. File: gsl-ref.info,  Node: Gamma Function,  Next: Gegenbauer Functions,  Prev: Fermi-Dirac Function,  Up: Special Functions
  201.  
  202. Gamma Function
  203. ==============
  204.  
  205.    /* Log[Gamma(x)], x not a negative integer  * Uses real Lanczos
  206. method.   * Returns the real part of Log[Gamma[x]] when x < 0,  * i.e.
  207. Log[|Gamma[x]|].   *  * exceptions: GSL_EDOM, GSL_EROUND  */ int
  208. gsl_sf_lngamma_impl(double x, gsl_sf_result * result); int
  209. gsl_sf_lngamma_e(double x, gsl_sf_result * result);
  210.  
  211.    /* Log[Gamma(x)], x not a negative integer  * Uses real Lanczos
  212. method.  Determines  * the sign of Gamma[x] as well as Log[|Gamma[x]|]
  213. for x < 0.   * So Gamma[x] = sgn * Exp[result_lg].   *  * exceptions:
  214. GSL_EDOM, GSL_EROUND  */ int gsl_sf_lngamma_sgn_impl(double x,
  215. gsl_sf_result * result_lg, double *sgn); int
  216. gsl_sf_lngamma_sgn_e(double x, gsl_sf_result * result_lg, double * sgn);
  217.  
  218.    /* Gamma(x), x not a negative integer  * Uses real Lanczos method.
  219. *  * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EROUND  */ int
  220. gsl_sf_gamma_impl(double x, gsl_sf_result * result); int
  221. gsl_sf_gamma_e(double x, gsl_sf_result * result);
  222.  
  223.    /* Regulated Gamma Function, x > 0  * Gamma^*(x) =
  224. Gamma(x)/(Sqrt[2Pi] x^(x-1/2) exp(-x))  *            = (1 + 1/(12x) +
  225. ...),  x->Inf  * A useful suggestion of Temme.   *  * exceptions:
  226. GSL_EDOM  */ int gsl_sf_gammastar_impl(double x, gsl_sf_result *
  227. result); int gsl_sf_gammastar_e(double x, gsl_sf_result * result);
  228.  
  229.    /* 1/Gamma(x)  * Uses real Lanczos method.   *  * exceptions:
  230. GSL_EUNDRFLW, GSL_EROUND  */ int gsl_sf_gammainv_impl(double x,
  231. gsl_sf_result * result); int gsl_sf_gammainv_e(double x, gsl_sf_result
  232. * result);
  233.  
  234.    /* Log[Gamma(z)] for z complex, z not a negative integer  * Uses
  235. complex Lanczos method.  Note that the phase part (arg)  * is not
  236. well-determined when |z| is very large, due  * to inevitable roundoff
  237. in restricting to (-Pi,Pi].   * This will raise the GSL_ELOSS exception
  238. when it occurs.   * The absolute value part (lnr), however, never
  239. suffers.   *  * Calculates:  *   lnr = log|Gamma(z)|  *   arg =
  240. arg(Gamma(z))  in (-Pi, Pi]  *  * exceptions: GSL_EDOM, GSL_ELOSS  */
  241. int gsl_sf_lngamma_complex_impl(double zr, double zi, gsl_sf_result *
  242. lnr, gsl_sf_result * arg); int gsl_sf_lngamma_complex_e(double zr,
  243. double zi, gsl_sf_result * lnr, gsl_sf_result * arg);
  244.  
  245.    /* x^n / n!   *  * x >= 0.0, n >= 0  * exceptions: GSL_EDOM,
  246. GSL_EOVRFLW, GSL_EUNDRFLW  */ int gsl_sf_taylorcoeff_impl(int n, double
  247. x, gsl_sf_result * result); int gsl_sf_taylorcoeff_e(int n, double x,
  248. gsl_sf_result * result);
  249.  
  250.    /* n!   *  * exceptions: GSL_EDOM, GSL_OVRFLW  */ int
  251. gsl_sf_fact_impl(unsigned int n, gsl_sf_result * result); int
  252. gsl_sf_fact_e(unsigned int n, gsl_sf_result * result);
  253.  
  254.    /* n!! = n(n-2)(n-4) ...   *  * exceptions: GSL_EDOM, GSL_OVRFLW  */
  255. int gsl_sf_doublefact_impl(unsigned int n, gsl_sf_result * result); int
  256. gsl_sf_doublefact_e(unsigned int n, gsl_sf_result * result);
  257.  
  258.    /* log(n!)   * Faster than ln(Gamma(n+1)) for n < 170; defers for
  259. larger n.   *  * exceptions: none  */ int gsl_sf_lnfact_impl(unsigned
  260. int n, gsl_sf_result * result); int gsl_sf_lnfact_e(unsigned int n,
  261. gsl_sf_result * result);
  262.  
  263.    /* log(n!!)   *  * exceptions: none  */ int
  264. gsl_sf_lndoublefact_impl(unsigned int n, gsl_sf_result * result); int
  265. gsl_sf_lndoublefact_e(unsigned int n, gsl_sf_result * result);
  266.  
  267.    /* log(n choose m)  *  * exceptions: GSL_EDOM  */ int
  268. gsl_sf_lnchoose_impl(unsigned int n, unsigned int m, gsl_sf_result *
  269. result); int gsl_sf_lnchoose_e(unsigned int n, unsigned int m,
  270. gsl_sf_result * result);
  271.  
  272.    /* n choose m  *  * exceptions: GSL_EDOM, GSL_EOVRFLW  */ int
  273. gsl_sf_choose_impl(unsigned int n, unsigned int m, gsl_sf_result *
  274. result); int gsl_sf_choose_e(unsigned int n, unsigned int m,
  275. gsl_sf_result * result);
  276.  
  277.    /* Logarithm of Pochammer (Apell) symbol  *   log( (a)_x )  *
  278. where (a)_x := Gamma[a + x]/Gamma[a]  *  * a > 0, a+x > 0  *  *
  279. exceptions:  GSL_EDOM  */ int gsl_sf_lnpoch_impl(double a, double x,
  280. gsl_sf_result * result); int gsl_sf_lnpoch_e(double a, double x,
  281. gsl_sf_result * result);
  282.  
  283.    /* Logarithm of Pochammer (Apell) symbol, with sign information.   *
  284.  result = log( |(a)_x| )  *   sgn    = sgn( (a)_x )  *   where (a)_x
  285. := Gamma[a + x]/Gamma[a]  *  * a != neg integer, a+x != neg integer  *
  286. * exceptions:  GSL_EDOM  */ int gsl_sf_lnpoch_sgn_impl(double a, double
  287. x, gsl_sf_result * result, double * sgn); int
  288. gsl_sf_lnpoch_sgn_e(double a, double x, gsl_sf_result * result, double
  289. * sgn);
  290.  
  291.    /* Pochammer (Apell) symbol  *   (a)_x := Gamma[a + x]/Gamma[x]  *
  292. * a != neg integer, a+x != neg integer  *  * exceptions:  GSL_EDOM,
  293. GSL_EOVRFLW  */ int gsl_sf_poch_impl(double a, double x, gsl_sf_result
  294. * result); int gsl_sf_poch_e(double a, double x, gsl_sf_result *
  295. result);
  296.  
  297.    /* Relative Pochammer (Apell) symbol  *   ((a,x) - 1)/x  *   where
  298. (a,x) = (a)_x := Gamma[a + x]/Gamma[a]  *  * exceptions:  GSL_EDOM  */
  299. int gsl_sf_pochrel_impl(double a, double x, gsl_sf_result * result);
  300. int gsl_sf_pochrel_e(double a, double x, gsl_sf_result * result);
  301.  
  302.    /* Normalized Incomplete Gamma Function  *  * Q(a,x) = 1/Gamma(a)
  303. Integral[ t^(a-1) e^(-t), (t,x,Infinity) ]  *  * a > 0, x >= 0  *  *
  304. exceptions: GSL_EDOM  */ int gsl_sf_gamma_inc_Q_impl(double a, double
  305. x, gsl_sf_result * result); int gsl_sf_gamma_inc_Q_e(double a, double
  306. x, gsl_sf_result * result);
  307.  
  308.    /* Complementary Normalized Incomplete Gamma Function  *  * P(a,x) =
  309. 1/Gamma(a) Integral[ t^(a-1) e^(-t), (t,0,x) ]  *  * a > 0, x >= 0  *
  310. * exceptions: GSL_EDOM  */ int gsl_sf_gamma_inc_P_impl(double a, double
  311. x, gsl_sf_result * result); int gsl_sf_gamma_inc_P_e(double a, double
  312. x, gsl_sf_result * result);
  313.  
  314.    /* Logarithm of Beta Function  * Log[B(a,b)]  *  * a > 0, b > 0  *
  315. exceptions: GSL_EDOM  */ int gsl_sf_lnbeta_impl(double a, double b,
  316. gsl_sf_result * result); int gsl_sf_lnbeta_e(double a, double b,
  317. gsl_sf_result * result);
  318.  
  319.    /* Beta Function  * B(a,b)  *  * a > 0, b > 0  * exceptions:
  320. GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW  */ int gsl_sf_beta_impl(double a,
  321. double b, gsl_sf_result * result); int gsl_sf_beta_e(double a, double
  322. b, gsl_sf_result * result);
  323.  
  324.    /* The maximum x such that gamma(x) is not  * considered an overflow.
  325. */ #define GSL_SF_GAMMA_XMAX  171.0
  326.  
  327. 
  328. File: gsl-ref.info,  Node: Gegenbauer Functions,  Next: Hypergeometric Functions,  Prev: Gamma Function,  Up: Special Functions
  329.  
  330. Gegenbauer Functions
  331. ====================
  332.  
  333.  - Function: int gsl_sf_gegenpoly_1_impl (double lambda, double x,
  334.           gsl_sf_result * result)
  335.  - Function: int gsl_sf_gegenpoly_2_impl (double lambda, double x,
  336.           gsl_sf_result * result)
  337.  - Function: int gsl_sf_gegenpoly_3_impl (double lambda, double x,
  338.           gsl_sf_result * result)
  339.  - Function: int gsl_sf_gegenpoly_1_e (double lambda, double x,
  340.           gsl_sf_result * result)
  341.  - Function: int gsl_sf_gegenpoly_2_e (double lambda, double x,
  342.           gsl_sf_result * result)
  343.  - Function: int gsl_sf_gegenpoly_3_e (double lambda, double x,
  344.           gsl_sf_result * result)
  345.      Evaluate Gegenbauer polynomials using explicit representations.
  346.      Exceptional Return Values: none
  347.  
  348.  - Function: int gsl_sf_gegenpoly_n_impl (int n, double lambda, double
  349.           x, gsl_sf_result * result)
  350.  - Function: int gsl_sf_gegenpoly_n_e (int n, double lambda, double x,
  351.           gsl_sf_result * result)
  352.      Evaluate Gegenbauer polynomials.  Domain: lambda > -1/2, n >= 0
  353.      Exceptional Return Values: GSL_EDOM
  354.  
  355.  - Function: int gsl_sf_gegenpoly_array_impl (int nmax, double lambda,
  356.           double x, double * result_array);
  357.  - Function: int gsl_sf_gegenpoly_array_e (int nmax, double lambda,
  358.           double x, double * result_array);
  359.      Calculate array of Gegenbauer polynomials.  Conditions: n = 0, 1,
  360.      2, ... nmax Domain: lambda > -1/2, nmax >= 0 Exceptional Return
  361.      Values: GSL_EDOM
  362.  
  363. 
  364. File: gsl-ref.info,  Node: Hypergeometric Functions,  Next: Laguerre Functions,  Prev: Gegenbauer Functions,  Up: Special Functions
  365.  
  366. Hypergeometric Functions
  367. ========================
  368.  
  369.    /* Hypergeometric function related to Bessel functions  * 0F1[c,x] =
  370. *            Gamma[c]    x^(1/2(1-c)) I_(c-1)(2 Sqrt[x])  *
  371. Gamma[c] (-x)^(1/2(1-c)) J_(c-1)(2 Sqrt[-x])  *  * exceptions:
  372. GSL_EOVRFLW, GSL_EUNDRFLW  */ int gsl_sf_hyperg_0F1_impl(double c,
  373. double x, gsl_sf_result * result); int gsl_sf_hyperg_0F1_e(double c,
  374. double x, gsl_sf_result * result);
  375.  
  376.    /* Confluent hypergeometric function  for integer parameters.   *
  377. 1F1[m,n,x] = M(m,n,x)  *  * exceptions:  */ int
  378. gsl_sf_hyperg_1F1_int_impl(int m, int n, double x, gsl_sf_result *
  379. result); int gsl_sf_hyperg_1F1_int_e(int m, int n, double x,
  380. gsl_sf_result * result);
  381.  
  382.    /* Confluent hypergeometric function.   * 1F1[a,b,x] = M(a,b,x)  *
  383. * exceptions:  */ int gsl_sf_hyperg_1F1_impl(double a, double b, double
  384. x, gsl_sf_result * result); int gsl_sf_hyperg_1F1_e(double a, double b,
  385. double x, gsl_sf_result * result);
  386.  
  387.    /* Confluent hypergeometric function for integer parameters.   *
  388. U(m,n,x)  *  * exceptions:  */ int gsl_sf_hyperg_U_int_impl(int m, int
  389. n, double x, gsl_sf_result * result); int gsl_sf_hyperg_U_int_e(int m,
  390. int n, double x, gsl_sf_result * result);
  391.  
  392.    /* Confluent hypergeometric function for integer parameters.   *
  393. U(m,n,x)  *  * exceptions:  */ int gsl_sf_hyperg_U_int_e10_impl(int m,
  394. int n, double x, gsl_sf_result_e10 * result); int
  395. gsl_sf_hyperg_U_int_e10_e(int m, int n, double x, gsl_sf_result_e10 *
  396. result);
  397.  
  398.    /* Confluent hypergeometric function.   * U(a,b,x)  *  * exceptions:
  399. */ int gsl_sf_hyperg_U_impl(double a, double b, double x,
  400. gsl_sf_result * result); int gsl_sf_hyperg_U_e(double a, double b,
  401. double x, gsl_sf_result * result);
  402.  
  403.    /* Confluent hypergeometric function.   * U(a,b,x)  *  * exceptions:
  404. */ int gsl_sf_hyperg_U_e10_impl(double a, double b, double x,
  405. gsl_sf_result_e10 * result); int gsl_sf_hyperg_U_e10_e(double a, double
  406. b, double x, gsl_sf_result_e10 * result);
  407.  
  408.    /* Gauss hypergeometric function 2F1[a,b,c,x]  * |x| < 1  *  *
  409. exceptions:  */ int gsl_sf_hyperg_2F1_impl(double a, double b, double
  410. c, double x, gsl_sf_result * result); int gsl_sf_hyperg_2F1_e(double a,
  411. double b, double c, double x, gsl_sf_result * result);
  412.  
  413.    /* Gauss hypergeometric function  * 2F1[aR + I aI, aR - I aI, c, x]
  414. * |x| < 1  *  * exceptions:  */ int gsl_sf_hyperg_2F1_conj_impl(double
  415. aR, double aI, double c, double x, gsl_sf_result * result); int
  416. gsl_sf_hyperg_2F1_conj_e(double aR, double aI, double c, double x,
  417. gsl_sf_result * result);
  418.  
  419.    /* Renormalized Gauss hypergeometric function  * 2F1[a,b,c,x] /
  420. Gamma[c]  * |x| < 1  *  * exceptions:  */ int
  421. gsl_sf_hyperg_2F1_renorm_impl(double a, double b, double c, double x,
  422. gsl_sf_result * result); int gsl_sf_hyperg_2F1_renorm_e(double a,
  423. double b, double c, double x, gsl_sf_result * result);
  424.  
  425.    /* Renormalized Gauss hypergeometric function  * 2F1[aR + I aI, aR -
  426. I aI, c, x] / Gamma[c]  * |x| < 1  *  * exceptions:  */ int
  427. gsl_sf_hyperg_2F1_conj_renorm_impl(double aR, double aI, double c,
  428. double x, gsl_sf_result * result); int
  429. gsl_sf_hyperg_2F1_conj_renorm_e(double aR, double aI, double c, double
  430. x, gsl_sf_result * result);
  431.  
  432.    /* Mysterious hypergeometric function.  The series representation  *
  433. is a divergent hypergeometric series.  However, for x < 0 we  * have
  434. 2F0(a,b,x) = (-1/x)^a U(a,1+a-b,-1/x)  *  * exceptions: GSL_EDOM  */
  435. int     gsl_sf_hyperg_2F0_impl(double a, double b, double x,
  436. gsl_sf_result * result); int     gsl_sf_hyperg_2F0_e(double a, double
  437. b, double x, gsl_sf_result * result);
  438.  
  439. 
  440. File: gsl-ref.info,  Node: Laguerre Functions,  Next: Legendre Functions and Spherical Harmonics,  Prev: Hypergeometric Functions,  Up: Special Functions
  441.  
  442. Laguerre Functions
  443. ==================
  444.  
  445.    The Laguerre polynomials are defined in terms of confluent
  446. hypergeometric functions as  L^a_n(x) = (a+1)_n / n! 1F1(-n,a+1,x) .
  447.  
  448.  - Function: int gsl_sf_laguerre_1_impl (double a, double x,
  449.           gsl_sf_result * result)
  450.  - Function: int gsl_sf_laguerre_2_impl (double a, double x,
  451.           gsl_sf_result * result)
  452.  - Function: int gsl_sf_laguerre_3_impl (double a, double x,
  453.           gsl_sf_result * result)
  454.  - Function: int gsl_sf_laguerre_1_e (double a, double x, gsl_sf_result
  455.           * result)
  456.  - Function: int gsl_sf_laguerre_2_e (double a, double x, gsl_sf_result
  457.           * result)
  458.  - Function: int gsl_sf_laguerre_3_e (double a, double x, gsl_sf_result
  459.           * result)
  460.      Evaluate generalized Laguerre polynomials using explicit
  461.      representations.  Exceptional Return Values: none
  462.  
  463.  - Function: int gsl_sf_laguerre_n_impl (const int n, const double a,
  464.           const double x, gsl_sf_result * result)
  465.  - Function: int gsl_sf_laguerre_n_e (int n, double a, double x,
  466.           gsl_sf_result * result)
  467.      Domain: a > -1.0, n >= 0 Evaluate generalized Laguerre polynomials.
  468.      Exceptional Return Values: GSL_EDOM
  469.  
  470. 
  471. File: gsl-ref.info,  Node: Legendre Functions and Spherical Harmonics,  Next: Logarithm and Related Functions,  Prev: Laguerre Functions,  Up: Special Functions
  472.  
  473. Legendre Functions and Spherical Harmonics
  474. ==========================================
  475.  
  476. Legendre Polynomials
  477. --------------------
  478.  
  479.    /* P_l(x), l=1,2,3  *
  480.  
  481.  - Function: int gsl_sf_legendre_P1_impl (double x, gsl_sf_result *
  482.           result)
  483.  - Function: int gsl_sf_legendre_P2_impl (double x, gsl_sf_result *
  484.           result)
  485.  - Function: int gsl_sf_legendre_P3_impl (double x, gsl_sf_result *
  486.           result)
  487.  - Function: int gsl_sf_legendre_P1_e (double x, gsl_sf_result * result)
  488.  - Function: int gsl_sf_legendre_P2_e (double x, gsl_sf_result * result)
  489.  - Function: int gsl_sf_legendre_P3_e (double x, gsl_sf_result * result)
  490.      Exceptional Return Values: none
  491.  
  492.    /* P_l(x)   l >= 0; |x| <= 1  *
  493.  
  494.  - Function: int gsl_sf_legendre_Pl_impl (int l, double x,
  495.           gsl_sf_result * result)
  496.  - Function: int gsl_sf_legendre_Pl_e (int l, double x, gsl_sf_result *
  497.           result)
  498.      Exceptional Return Values: GSL_EDOM
  499.  
  500.    /* P_l(x) for l=0,...,lmax; |x| <= 1  *
  501.  
  502.  - Function: int gsl_sf_legendre_Pl_array_impl (int lmax, double x,
  503.           double * result_array)
  504.  - Function: int gsl_sf_legendre_Pl_array_e (int lmax, double x, double
  505.           * result_array)
  506.      Exceptional Return Values: GSL_EDOM
  507.  
  508.    /* Q_0(x), x > -1, x != 1  *
  509.  
  510.  - Function: int gsl_sf_legendre_Q0_impl (double x, gsl_sf_result *
  511.           result)
  512.  - Function: int gsl_sf_legendre_Q0_e (double x, gsl_sf_result * result)
  513.      Exceptional Return Values: GSL_EDOM
  514.  
  515.    /* Q_1(x), x > -1, x != 1  *
  516.  
  517.  - Function: int gsl_sf_legendre_Q1_impl (double x, gsl_sf_result *
  518.           result)
  519.  - Function: int gsl_sf_legendre_Q1_e (double x, gsl_sf_result * result)
  520.      Exceptional Return Values: GSL_EDOM
  521.  
  522.    /* Q_l(x), x > -1, x != 1, l >= 0  *
  523.  
  524.  - Function: int gsl_sf_legendre_Ql_impl (int l, double x,
  525.           gsl_sf_result * result)
  526.  - Function: int gsl_sf_legendre_Ql_e (int l, double x, gsl_sf_result *
  527.           result)
  528.      Exceptional Return Values: GSL_EDOM
  529.  
  530. Associated Legendre Polynomials and Spherical Harmonics
  531. -------------------------------------------------------
  532.  
  533.    /* P_l^m(x)  m >= 0; l >= m; |x| <= 1.0  *  * Note that this
  534. function grows combinatorially with l.   * Therefore we can easily
  535. generate an overflow for l larger  * than about 150.   *  * There is no
  536. trouble for small m, but when m and l are both large,  * then there
  537. will be trouble.  Rather than allow overflows, these  * functions
  538. refuse to calculate when they can sense that l and m are  * too big.   *
  539. * If you really want to calculate a spherical harmonic, then DO NOT  *
  540. use this.  Instead use legendre_sphPlm() below, which  uses a similar
  541. * recursion, but with the normalized functions.   *
  542.  
  543.  - Function: int gsl_sf_legendre_Plm_impl (int l, int m, double x,
  544.           gsl_sf_result * result)
  545.  - Function: int gsl_sf_legendre_Plm_e (int l, int m, double x,
  546.           gsl_sf_result * result)
  547.      Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW
  548.  
  549.    /* P_l^m(x)  m >= 0; l >= m; |x| <= 1.0  * l=|m|,...,lmax  *
  550.  
  551.  - Function: int gsl_sf_legendre_Plm_array_impl (int lmax, int m,
  552.           double x, double * result_array)
  553.  - Function: int gsl_sf_legendre_Plm_array_e (int lmax, int m, double
  554.           x, double * result_array)
  555.      Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW
  556.  
  557.    /* P_l^m(x), normalized properly for use in spherical harmonics  * m
  558. >= 0; l >= m; |x| <= 1.0  *  * There is no overflow problem, as there
  559. is for the  * standard normalization of P_l^m(x).   *  * Specifically,
  560. it returns:  *  *        sqrt((2l+1)/(4pi)) sqrt((l-m)!/(l+m)!) P_l^m(x)
  561. *
  562.  
  563.  - Function: int gsl_sf_legendre_sphPlm_impl (int l, int m, double x,
  564.           gsl_sf_result * result)
  565.  - Function: int gsl_sf_legendre_sphPlm_e (int l, int m, double x,
  566.           gsl_sf_result * result)
  567.      Exceptional Return Values: GSL_EDOM
  568.  
  569.    /* P_l^m(x), normalized properly for use in spherical harmonics  * m
  570. >= 0; l >= m; |x| <= 1.0  * l=|m|,...,lmax  *
  571.  
  572.  - Function: int gsl_sf_legendre_sphPlm_array_impl (int lmax, int m,
  573.           double x, double * result_array)
  574.  - Function: int gsl_sf_legendre_sphPlm_array_e (int lmax, int m,
  575.           double x, double * result_array)
  576.      Exceptional Return Values: GSL_EDOM
  577.  
  578.    /* size of result_array[] needed for the array versions of Plm  *
  579. (lmax - m + 1)  */
  580.  
  581.  - Function: int gsl_sf_legendre_array_size (const int lmax, const int
  582.           m)
  583.      Exceptional Return Values: none
  584.  
  585. Conical Functions
  586. -----------------
  587.  
  588.    /* Irregular Spherical Conical Function  * P^(1/2)_(-1/2 + I
  589. lambda)(x)  *  * x > -1.0
  590.  
  591.  - Function: int gsl_sf_conicalP_half_impl (double lambda, double x,
  592.           gsl_sf_result * result)
  593.  - Function: int gsl_sf_conicalP_half_e (double lambda, double x,
  594.           gsl_sf_result * result)
  595.      Exceptional Return Values: GSL_EDOM
  596.  
  597.    /* Regular Spherical Conical Function  * P^(-1/2)_(-1/2 + I
  598. lambda)(x)  *  * x > -1.0
  599.  
  600.  - Function: int gsl_sf_conicalP_mhalf_impl (double lambda, double x,
  601.           gsl_sf_result * result)
  602.  - Function: int gsl_sf_conicalP_mhalf_e (double lambda, double x,
  603.           gsl_sf_result * result)
  604.      Exceptional Return Values: GSL_EDOM
  605.  
  606.    /* Conical Function  * P^(0)_(-1/2 + I lambda)(x)  *  * x > -1.0
  607.  
  608.  - Function: int gsl_sf_conicalP_0_impl (double lambda, double x,
  609.           gsl_sf_result * result)
  610.  - Function: int gsl_sf_conicalP_0_e (double lambda, double x,
  611.           gsl_sf_result * result)
  612.      Exceptional Return Values: GSL_EDOM
  613.  
  614.    /* Conical Function  * P^(1)_(-1/2 + I lambda)(x)  *  * x > -1.0
  615.  
  616.  - Function: int gsl_sf_conicalP_1_impl (double lambda, double x,
  617.           gsl_sf_result * result)
  618.  - Function: int gsl_sf_conicalP_1_e (double lambda, double x,
  619.           gsl_sf_result * result)
  620.      Exceptional Return Values: GSL_EDOM
  621.  
  622.    /* Regular Spherical Conical Function  * P^(-1/2-l)_(-1/2 + I
  623. lambda)(x)  *  * x > -1.0, l >= -1
  624.  
  625.  - Function: int gsl_sf_conicalP_sph_reg_impl (int l, double lambda,
  626.           double x, gsl_sf_result * result)
  627.  - Function: int gsl_sf_conicalP_sph_reg_e (int l, double lambda,
  628.           double x, gsl_sf_result * result)
  629.      Exceptional Return Values: GSL_EDOM
  630.  
  631.    /* Regular Cylindrical Conical Function  * P^(-m)_(-1/2 + I
  632. lambda)(x)  *  * x > -1.0, m >= -1
  633.  
  634.  - Function: int gsl_sf_conicalP_cyl_reg_impl (int m, double lambda,
  635.           double x, gsl_sf_result * result)
  636.  - Function: int gsl_sf_conicalP_cyl_reg_e (int m, double lambda,
  637.           double x, gsl_sf_result * result)
  638.      Exceptional Return Values: GSL_EDOM
  639.  
  640. Radial Functions for Hyperbolic Space
  641. -------------------------------------
  642.  
  643.    /* The following spherical functions are specializations  * of
  644. Legendre functions which give the regular eigenfunctions  * of the
  645. Laplacian on a 3-dimensional hyperbolic space.   * Of particular
  646. interest is the flat limit, which is  * Flat-Lim := (lambda->Inf,
  647. eta->0, lambda*eta fixed).   */
  648.  
  649.    /* Zeroth radial eigenfunction of the Laplacian on the  *
  650. 3-dimensional hyperbolic space.   *  * legendre_H3d_0(lambda,eta) :=
  651. sin(lambda*eta)/(lambda*sinh(eta))  *  * Normalization:  * Flat-Lim
  652. legendre_H3d_0(lambda,eta) = j_0(lambda*eta)  *  * eta >= 0.0
  653.  
  654.  - Function: int gsl_sf_legendre_H3d_0_impl (double lambda, double eta,
  655.           gsl_sf_result * result)
  656.  - Function: int gsl_sf_legendre_H3d_0_e (double lambda, double eta,
  657.           gsl_sf_result * result)
  658.      Exceptional Return Values: GSL_EDOM
  659.  
  660.    /* First radial eigenfunction of the Laplacian on the  *
  661. 3-dimensional hyperbolic space.   *  * legendre_H3d_1(lambda,eta) :=  *
  662.   1/sqrt(lambda^2 + 1) sin(lam eta)/(lam sinh(eta))  *    (coth(eta) -
  663. lambda cot(lambda*eta))  *  * Normalization:  * Flat-Lim
  664. legendre_H3d_1(lambda,eta) = j_1(lambda*eta)  *  * eta >= 0.0
  665.  
  666.  - Function: int gsl_sf_legendre_H3d_1_impl (double lambda, double eta,
  667.           gsl_sf_result * result)
  668.  - Function: int gsl_sf_legendre_H3d_1_e (double lambda, double eta,
  669.           gsl_sf_result * result)
  670.      Exceptional Return Values: GSL_EDOM
  671.  
  672.    /* l'th radial eigenfunction of the Laplacian on the  *
  673. 3-dimensional hyperbolic space.   *  * Normalization:  * Flat-Lim
  674. legendre_H3d_l(l,lambda,eta) = j_l(lambda*eta)  *  * eta >= 0.0, l >= 0
  675.  
  676.  - Function: int gsl_sf_legendre_H3d_impl (int l, double lambda, double
  677.           eta, gsl_sf_result * result)
  678.  - Function: int gsl_sf_legendre_H3d_e (int l, double lambda, double
  679.           eta, gsl_sf_result * result)
  680.      Exceptional Return Values: GSL_EDOM
  681.  
  682.    /* Array of H3d(ell),  0 <= ell <= lmax  */
  683.  
  684.  - Function: int gsl_sf_legendre_H3d_array_impl (int lmax, double
  685.           lambda, double eta, double * result_array)
  686.  - Function: int gsl_sf_legendre_H3d_array_e (int lmax, double lambda,
  687.           double eta, double * result_array)
  688.      Exceptional Return Values:
  689.  
  690. 
  691. File: gsl-ref.info,  Node: Logarithm and Related Functions,  Next: Polynomial Manipulation,  Prev: Legendre Functions and Spherical Harmonics,  Up: Special Functions
  692.  
  693. Logarithm and Related Functions
  694. ===============================
  695.  
  696.  - Function: int gsl_sf_log_impl (double x, gsl_sf_result * result)
  697.  - Function: int gsl_sf_log_e (double x, gsl_sf_result * result)
  698.      Domain: x > 0.0 Exceptional Return Values: GSL_EDOM
  699.  
  700.  - Function: int gsl_sf_log_abs_impl (double x, gsl_sf_result * result)
  701.  - Function: int gsl_sf_log_abs_e (double x, gsl_sf_result * result)
  702.      \log(|x|) Domain: x != 0.0 Exceptional Return Values: GSL_EDOM
  703.  
  704.    /* Complex Logarithm  *   exp(lnr + I theta) = zr + I zi  * Returns
  705. argument in [-pi,pi].   *
  706.  
  707.  - Function: int gsl_sf_complex_log_impl (double zr, double zi,
  708.           gsl_sf_result * lnr, gsl_sf_result * theta)
  709.  - Function: int gsl_sf_complex_log_e (double zr, double zi,
  710.           gsl_sf_result * lnr, gsl_sf_result * theta)
  711.      Exceptional Return Values: GSL_EDOM
  712.  
  713.  - Function: int gsl_sf_log_1plusx_impl (double x, gsl_sf_result *
  714.           result)
  715.  - Function: int gsl_sf_log_1plusx_e (double x, gsl_sf_result * result)
  716.      \log(1 + x) Domain: x > -1.0 Exceptional Return Values: GSL_EDOM
  717.  
  718.  - Function: int gsl_sf_log_1plusx_mx_impl (double x, gsl_sf_result *
  719.           result)
  720.  - Function: int gsl_sf_log_1plusx_mx_e (double x, gsl_sf_result *
  721.           result)
  722.      \log(1 + x) - x Domain: x > -1.0 Exceptional Return Values:
  723.      GSL_EDOM
  724.  
  725. 
  726. File: gsl-ref.info,  Node: Polynomial Manipulation,  Next: Power Function,  Prev: Logarithm and Related Functions,  Up: Special Functions
  727.  
  728. Polynomial Manipulation
  729. =======================
  730.  
  731.  - Function: double gsl_sf_poly_eval (const double c[], const int len,
  732.           const double x)
  733.      Evaluate a polynomial,  c[0] + c[1] x + c[2] x^2 + ... + c[len-1]
  734.      x^(len-1) .  Exceptional Return Values: none
  735.  
  736. 
  737. File: gsl-ref.info,  Node: Power Function,  Next: Psi (Digamma) Function,  Prev: Polynomial Manipulation,  Up: Special Functions
  738.  
  739. Power Function
  740. ==============
  741.  
  742.    A common complaint about the standard C library is its lack of a
  743. function for calculating (small) integer powers. GSL provides a simple
  744. function to fill this gap.  For reasons of efficiency, these functions
  745. do not check for overflow or underflow conditions.
  746.  
  747.    int     gsl_sf_pow_int_impl(double x, int n, gsl_sf_result * result);
  748. int     gsl_sf_pow_int_e(double x, int n, gsl_sf_result * result);
  749.  
  750.  - Function: double gsl_sf_pow_int (double x, int n)
  751.  
  752.    double gsl_sf_pow_2(const double x); double gsl_sf_pow_3(const
  753. double x); double gsl_sf_pow_4(const double x); double
  754. gsl_sf_pow_5(const double x); double gsl_sf_pow_6(const double x);
  755. double gsl_sf_pow_7(const double x); double gsl_sf_pow_8(const double
  756. x); double gsl_sf_pow_9(const double x);
  757.  
  758.      #include <gsl/gsl_sf_pow_int.h>
  759.      double y = gsl_sf_pow_int(3.0, 12)
  760.  
  761. 
  762. File: gsl-ref.info,  Node: Psi (Digamma) Function,  Next: Synchrotron Functions,  Prev: Power Function,  Up: Special Functions
  763.  
  764. Psi (Digamma) Function
  765. ======================
  766.  
  767.    The poly-gamma functions are defined by  psi(m,x) = (d/dx)^m
  768. psi(0,x) = (d/dx)^(m+1) log(gamma(x)) , where  \psi(0,x)  is called the
  769. di-gamma function.
  770.  
  771. Digamma Function
  772. ----------------
  773.  
  774.  - Function: int gsl_sf_psi_int_impl (int n, gsl_sf_result * result)
  775.  - Function: int gsl_sf_psi_int_e (int n, gsl_sf_result * result)
  776.      Domain: n integer, n > 0 Exceptional Return Values: GSL_EDOM
  777.  
  778.  - Function: int gsl_sf_psi_impl (double x, gsl_sf_result * result)
  779.  - Function: int gsl_sf_psi_e (double x, gsl_sf_result * result)
  780.      Domain: x != 0.0 Exceptional Return Values: GSL_EDOM, GSL_ELOSS
  781.  
  782.    /* Di-Gamma Function Re[psi(1 + I y)]  *  * exceptions: none  */
  783.  
  784.  - Function: int gsl_sf_psi_1piy_impl (double y, gsl_sf_result * result)
  785.  - Function: int gsl_sf_psi_1piy_e (double y, gsl_sf_result * result)
  786.      Computes  Re psi(1 + I y) .  Exceptional Return Values: none
  787.  
  788. Trigamma Function
  789. -----------------
  790.  
  791.  - Function: int gsl_sf_psi_1_int_impl (int n, gsl_sf_result * result)
  792.  - Function: int gsl_sf_psi_1_int_e (int n, gsl_sf_result * result)
  793.      Domain: n integer, n > 0 Exceptional Return Values: GSL_EDOM
  794.  
  795. Polygamma Function
  796. ------------------
  797.  
  798.  - Function: int gsl_sf_psi_n_impl (int m, double x, gsl_sf_result *
  799.           result)
  800.  - Function: int gsl_sf_psi_n_e (int m, double x, gsl_sf_result *
  801.           result)
  802.      Domain: m >= 0, x > 0.0 Exceptional Return Values: GSL_EDOM
  803.  
  804. 
  805. File: gsl-ref.info,  Node: Synchrotron Functions,  Next: Transport Functions,  Prev: Psi (Digamma) Function,  Up: Special Functions
  806.  
  807. Synchrotron Functions
  808. =====================
  809.  
  810.    The first synchrotron function is defined by the integral
  811. representation  synchrotron_1(x) = x \int_x^\infty dt K_(5/3)(t) .
  812.  
  813.  - Function: int gsl_sf_synchrotron_1_impl (double x, gsl_sf_result *
  814.           result)
  815.  - Function: int gsl_sf_synchrotron_1_e (double x, gsl_sf_result *
  816.           result)
  817.      Domain: x >= 0.0 Exceptional Return Values: GSL_EDOM, GSL_EUNDRFLW
  818.  
  819.    The second synchroton function is defined by  synchrotron_2(x) = x *
  820. K_(2/3)(x) .
  821.  
  822.  - Function: int gsl_sf_synchrotron_2_impl (double x, gsl_sf_result *
  823.           result)
  824.  - Function: int gsl_sf_synchrotron_2_e (double x, gsl_sf_result *
  825.           result)
  826.      Domain: x >= 0.0 Exceptional Return Values: GSL_EDOM, GSL_EUNDRFLW
  827.  
  828. 
  829. File: gsl-ref.info,  Node: Transport Functions,  Next: Trigonometric Functions,  Prev: Synchrotron Functions,  Up: Special Functions
  830.  
  831. Transport Functions
  832. ===================
  833.  
  834.    The transport functions are defined by the integral representations
  835. J(n,x) := \int_0^x dt t^n e^t /(e^t - 1)^2 .
  836.  
  837.  - Function: int gsl_sf_transport_2_impl (double x, gsl_sf_result *
  838.           result)
  839.  - Function: int gsl_sf_transport_2_e (double x, gsl_sf_result * result)
  840.      Calculates J(2,x).  Exceptional Return Values: GSL_EDOM
  841.  
  842.  - Function: int gsl_sf_transport_3_impl (double x, gsl_sf_result *
  843.           result)
  844.  - Function: int gsl_sf_transport_3_e (double x, gsl_sf_result * result)
  845.      Calculates J(3,x).  Exceptional Return Values: GSL_EDOM,
  846.      GSL_EUNDRFLW
  847.  
  848.  - Function: int gsl_sf_transport_4_impl (double x, gsl_sf_result *
  849.           result)
  850.  - Function: int gsl_sf_transport_4_e (double x, gsl_sf_result * result)
  851.      Calculates J(4,x).  Exceptional Return Values: GSL_EDOM,
  852.      GSL_EUNDRFLW
  853.  
  854.  - Function: int gsl_sf_transport_5_impl (double x, gsl_sf_result *
  855.           result)
  856.  - Function: int gsl_sf_transport_5_e (double x, gsl_sf_result * result)
  857.      Calculates J(5,x).  Exceptional Return Values: GSL_EDOM,
  858.      GSL_EUNDRFLW
  859.  
  860. 
  861. File: gsl-ref.info,  Node: Trigonometric Functions,  Next: Zeta Functions,  Prev: Transport Functions,  Up: Special Functions
  862.  
  863. Trigonometric Functions
  864. =======================
  865.  
  866. Trigonometric Functions
  867. -----------------------
  868.  
  869.    These simple trogonometric functions are important because we want
  870. to control the error estimate, and trying to guess the error for the
  871. standard library implementation every time it is used would be a little
  872. goofy.
  873.  
  874.  - Function: int gsl_sf_sin_impl (double x, gsl_sf_result * result);
  875.  - Function: int gsl_sf_sin_e (double x, gsl_sf_result * result);
  876.      Exceptional Return Values:
  877.  
  878.  - Function: int gsl_sf_cos_impl (double x, gsl_sf_result * result);
  879.  - Function: int gsl_sf_cos_e (double x, gsl_sf_result * result);
  880.      Exceptional Return Values:
  881.  
  882.  - Function: int gsl_sf_hypot_impl (double x, double y, gsl_sf_result *
  883.           result);
  884.  - Function: int gsl_sf_hypot_e (double x, double y, gsl_sf_result *
  885.           result);
  886.      Exceptional Return Values:
  887.  
  888.  - Function: int gsl_sf_complex_sin_impl (double zr, double zi,
  889.           gsl_sf_result * szr, gsl_sf_result * szi);
  890.  - Function: int gsl_sf_complex_sin_e (double zr, double zi,
  891.           gsl_sf_result * szr, gsl_sf_result * szi);
  892.      Sin(z) Exceptional Return Values: GSL_EOVRFLW
  893.  
  894.  - Function: int gsl_sf_complex_cos_impl (double zr, double zi,
  895.           gsl_sf_result * czr, gsl_sf_result * czi);
  896.  - Function: int gsl_sf_complex_cos_e (double zr, double zi,
  897.           gsl_sf_result * czr, gsl_sf_result * czi);
  898.      Cos(z) Exceptional Return Values: GSL_EOVRFLW
  899.  
  900.  - Function: int gsl_sf_complex_logsin_impl (double zr, double zi,
  901.           gsl_sf_result * lszr, gsl_sf_result * lszi);
  902.  - Function: int gsl_sf_complex_logsin_e (double zr, double zi,
  903.           gsl_sf_result * lszr, gsl_sf_result * lszi);
  904.      Log(Sin(z)) Exceptional Return Values: GSL_EDOM, GSL_ELOSS
  905.  
  906.  - Function: int gsl_sf_sinc_impl (double x, gsl_sf_result * result);
  907.  - Function: int gsl_sf_sinc_e (double x, gsl_sf_result * result);
  908.      Sinc(x) = sin(pi x) / (pi x) Exceptional Return Values: none
  909.  
  910.  - Function: int gsl_sf_lnsinh_impl (double x, gsl_sf_result * result);
  911.  - Function: int gsl_sf_lnsinh_e (double x, gsl_sf_result * result);
  912.      Log(Sinh(x)) Domain: x > 0 Exceptional Return Values: GSL_EDOM
  913.  
  914.  - Function: int gsl_sf_lncosh_impl (double x, gsl_sf_result * result);
  915.  - Function: int gsl_sf_lncosh_e (double x, gsl_sf_result * result);
  916.      Log(Cosh(x)) Exceptional Return Values: none
  917.  
  918. Conversion Functions
  919. --------------------
  920.  
  921.  - Function: int gsl_sf_polar_to_rect_impl (double r, double theta,
  922.           gsl_sf_result * x, gsl_sf_result * y);
  923.  - Function: int gsl_sf_polar_to_rect_e (double r, double theta,
  924.           gsl_sf_result * x, gsl_sf_result * y);
  925.      Convert polar to rectlinear coordinates.  Exceptional Return
  926.      Values: GSL_ELOSS
  927.  
  928.  - Function: int gsl_sf_rect_to_polar_impl (double x, double y,
  929.           gsl_sf_result * r, gsl_sf_result * theta)
  930.  - Function: int gsl_sf_rect_to_polar_e (double x, double y,
  931.           gsl_sf_result * r, gsl_sf_result * theta)
  932.      Convert rectilinear to polar coordinates.  Return argument in
  933.      range [-pi, pi].  Exceptional Return Values: GSL_EDOM
  934.  
  935. Restriction Functions
  936. ---------------------
  937.  
  938.  - Function: int gsl_sf_angle_restrict_symm_impl (double * theta);
  939.  - Function: int gsl_sf_angle_restrict_symm_e (double * theta);
  940.      Force an angle to lie in the range (-pi,pi].  Exceptional Return
  941.      Values: GSL_ELOSS
  942.  
  943.  - Function: int gsl_sf_angle_restrict_pos_impl (double * theta);
  944.  - Function: int gsl_sf_angle_restrict_pos_e (double * theta);
  945.      Force an angle to lie in the range [0, 2pi).  Exceptional Return
  946.      Values: GSL_ELOSS
  947.  
  948.    Trigonometric Functions With Error Estimate
  949.  
  950.  - Function: int gsl_sf_sin_err_impl (double x, double dx,
  951.           gsl_sf_result * result);
  952.  - Function: int gsl_sf_sin_err_e (double x, double dx, gsl_sf_result *
  953.           result);
  954.  
  955.  - Function: int gsl_sf_cos_err_impl (double x, double dx,
  956.           gsl_sf_result * result);
  957.  - Function: int gsl_sf_cos_err_e (double x, double dx, gsl_sf_result *
  958.           result);
  959.  
  960. 
  961. File: gsl-ref.info,  Node: Zeta Functions,  Prev: Trigonometric Functions,  Up: Special Functions
  962.  
  963. Zeta Functions
  964. ==============
  965.  
  966. Riemann Zeta Function
  967. ---------------------
  968.  
  969.    The Riemann zeta function is defined by
  970.  
  971.  - Function: int gsl_sf_zeta_int_impl (int n, gsl_sf_result * result)
  972.  - Function: int gsl_sf_zeta_int_e (int n, gsl_sf_result * result)
  973.      Domain: n integer, n != 1 Exceptional Return Values: GSL_EDOM,
  974.      GSL_EOVRFLW
  975.  
  976.  - Function: int gsl_sf_zeta_impl (double s, gsl_sf_result * result)
  977.  - Function: int gsl_sf_zeta_e (double s, gsl_sf_result * result)
  978.      Domain: s != 1.0 Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW
  979.  
  980. Hurwitz Zeta Function
  981. ---------------------
  982.  
  983.    The Hurwitz zeta function is defined by  zeta(s,q) = \sum_0^\infty
  984. (k+q)^(-s) .
  985.  
  986.  - Function: int gsl_sf_hzeta_impl (double s, double q, gsl_sf_result *
  987.           result)
  988.  - Function: int gsl_sf_hzeta_e (double s, double q, gsl_sf_result *
  989.           result)
  990.      Domain: s > 1.0, q > 0.0 Exceptional Return Values: GSL_EDOM,
  991.      GSL_EUNDRFLW, GSL_EOVRFLW
  992.  
  993. Eta Function
  994. ------------
  995.  
  996.    The eta function is defined by  \eta(s) = (1-2^(1-s)) zeta(s) .
  997.  
  998.  - Function: int gsl_sf_eta_int_impl (int n, gsl_sf_result * result)
  999.  - Function: int gsl_sf_eta_int_e (int n, gsl_sf_result * result)
  1000.      Exceptional Return Values: GSL_EUNDRFLW, GSL_EOVRFLW
  1001.  
  1002.  - Function: int gsl_sf_eta_impl (double s, gsl_sf_result * result)
  1003.  - Function: int gsl_sf_eta_e (double s, gsl_sf_result * result)
  1004.      Exceptional Return Values: GSL_EUNDRFLW, GSL_EOVRFLW
  1005.  
  1006. 
  1007. File: gsl-ref.info,  Node: Roots of Polynomials,  Next: Interpolation,  Prev: Special Functions,  Up: Top
  1008.  
  1009. Roots of Polynomials
  1010. ********************
  1011.  
  1012.    This chapter describes functions for solving polynomial equations.
  1013. There are routines for finding real and complex roots of quadratic and
  1014. cubic equations using analytic methods.  An iterative polynomial solver
  1015. is also available for finding the roots of general real polynomials (of
  1016. any order).  The functions are defined in the header file `gsl_poly.h'.
  1017.  
  1018. * Menu:
  1019.  
  1020. * Quadratic equations::
  1021. * Cubic equations::
  1022. * General polynomial equations::
  1023. * Roots of Polynomials Examples::
  1024. * Roots of Polynomials References and Further Reading::
  1025.  
  1026. 
  1027. File: gsl-ref.info,  Node: Quadratic equations,  Next: Cubic equations,  Up: Roots of Polynomials
  1028.  
  1029. Quadratic equations
  1030. ===================
  1031.  
  1032.  - Function: int gsl_poly_solve_quadratic (double A, double B, double
  1033.           C, double *X0, double *X1)
  1034.      This function finds the real roots of the quadratic equation,
  1035.  
  1036.           a x^2 + b x + c = 0
  1037.  
  1038.      The number of real roots (either zero or two) is returned, and
  1039.      their locations are stored in X0 and X1.  If no real roots are
  1040.      found then X0 and X1 are not modified.  When two real roots are
  1041.      found they are stored in X0 and X1 in ascending order.  The case
  1042.      of coincident roots is not considered special.  For example
  1043.      (x-1)^2=0 will have two roots, which happen to have exactly equal
  1044.      values.
  1045.  
  1046.      The number of roots found depends on the sign of the discriminant
  1047.      b^2 - 4 a c.  This will be subject to rounding and cancellation
  1048.      errors when computed in double precision, and will also be subject
  1049.      to errors if the coefficients of the polynomial are inexact.
  1050.      These errors may cause a discrete change in the number of roots.
  1051.      However, for polynomials with small integer coefficients the
  1052.      discriminant can always be computed exactly.
  1053.  
  1054.  
  1055.  - Function: int gsl_poly_complex_solve_quadratic (double A, double B,
  1056.           double C, gsl_complex *Z0, gsl_complex *Z1)
  1057.      This function finds the complex roots of the quadratic equation,
  1058.  
  1059.           a z^2 + b z + c = 0
  1060.  
  1061.      The number of complex roots is returned (always two) and the
  1062.      locations of the roots are stored in Z0 and Z1.  The roots are
  1063.      returned in ascending order, sorted first by their real components
  1064.      and then by their imaginary components.
  1065.  
  1066.  
  1067. 
  1068. File: gsl-ref.info,  Node: Cubic equations,  Next: General polynomial equations,  Prev: Quadratic equations,  Up: Roots of Polynomials
  1069.  
  1070. Cubic equations
  1071. ===============
  1072.  
  1073.  - Function: int gsl_poly_solve_cubic (double A, double B, double C,
  1074.           double *X0, double *X1, double *X2)
  1075.      This function finds the real roots of the cubic equation,
  1076.  
  1077.           x^3 + a x^2 + b x + c = 0
  1078.  
  1079.      with a leading coefficient of unity.  The number of real roots
  1080.      (either one or three) is returned, and their locations are stored
  1081.      in X0, X1 and X2.  If one real root is found then only X0 is
  1082.      modified.  When three real roots are found they are stored in X0,
  1083.      X1 and X2 in ascending order.  The case of coincident roots is not
  1084.      considered special.  For example, the equation (x-1)^3=0 will have
  1085.      three roots with exactly equal values.
  1086.  
  1087.  
  1088.  - Function: int gsl_poly_complex_solve_cubic (double A, double B,
  1089.           double C, gsl_complex *Z0, gsl_complex *Z1, gsl_complex *Z2)
  1090.      This function finds the complex roots of the cubic equation,
  1091.  
  1092.           z^3 + a z^2 + b z + c = 0
  1093.  
  1094.      The number of complex roots is returned (always three) and the
  1095.      locations of the roots are stored in Z0, Z1 and Z2.  The roots are
  1096.      returned in ascending order, sorted first by their real components
  1097.      and then by their imaginary components.
  1098.  
  1099.  
  1100. 
  1101. File: gsl-ref.info,  Node: General polynomial equations,  Next: Roots of Polynomials Examples,  Prev: Cubic equations,  Up: Roots of Polynomials
  1102.  
  1103. General polynomial equations
  1104. ============================
  1105.  
  1106.    The roots of polynomial equations cannot be found analytically beyond
  1107. the special cases of the quadratic, cubic and quartic equation.  The
  1108. algorithm described in this section uses an iterative method to find the
  1109. approximate locations of roots of higher order polynomials.
  1110.  
  1111.  - Function: gsl_poly_complex_workspace *
  1112. gsl_poly_complex_workspace_alloc (size_t N)
  1113.      This function allocates space for a `gsl_poly_complex_workspace'
  1114.      struct and a scratch area suitable for solving a polynomial with N
  1115.      coefficients using the routine `gsl_poly_complex_solve'.
  1116.  
  1117.      The function returns a pointer to the newly allocated
  1118.      `gsl_poly_complex_workspace' if no errors were detected, and 0 in
  1119.      the case of error.
  1120.  
  1121.  - Function: void gsl_poly_complex_workspace_free
  1122.           (gsl_poly_complex_workspace * W)
  1123.      This function frees all the memory associated with the workspace W.
  1124.  
  1125.  - Function: int gsl_poly_complex_solve (const double * A, size_t N,
  1126.           gsl_poly_complex_workspace * W, gsl_complex_packed_ptr Z)
  1127.      This function computes the roots of the general polynomial P(x) =
  1128.      a_0 + a_1 x + a_2 x^2 + ... + a_{n-1} x^{n-1} using balanced-QR
  1129.      reduction of the companion matrix.  The parameter N specifies the
  1130.      length of the coefficient array.  The coefficient of the highest
  1131.      order term must be non-zero.  The function requires a workspace W
  1132.      of the appropriate size.  The n-1 roots are returned in the packed
  1133.      complex array Z of length 2(n-1), alternating real and imaginary
  1134.      parts.
  1135.  
  1136.      The functions returns `GSL_SUCCESS' if all the roots are found and
  1137.      `GSL_EFAILED' if the QR reduction does not converge.
  1138.  
  1139. 
  1140. File: gsl-ref.info,  Node: Roots of Polynomials Examples,  Next: Roots of Polynomials References and Further Reading,  Prev: General polynomial equations,  Up: Roots of Polynomials
  1141.  
  1142. Examples
  1143. ========
  1144.  
  1145.    To demonstrate the use of the general polynomial solver we will take
  1146. the polynomial P(x) = x^5 - 1 which has the following roots,
  1147.  
  1148.      1, e^{2\pi i /5}, e^{4\pi i /5}, e^{6\pi i /5}, e^{8\pi i /5}
  1149.  
  1150. The following program will find these roots.
  1151.  
  1152.      #include <stdio.h>
  1153.      #include <gsl/gsl_poly.h>
  1154.      
  1155.      int
  1156.      main ()
  1157.      {
  1158.        int i;
  1159.        double a[6] = { -1, 0, 0, 0, 0, 1 };  /*  P(x) =  -1 + x^5  */
  1160.        double z[10];
  1161.      
  1162.        gsl_poly_complex_workspace * w = gsl_poly_complex_workspace_alloc (6) ;
  1163.      
  1164.        gsl_poly_complex_solve (a, 6, w, z) ;
  1165.      
  1166.        gsl_poly_complex_workspace_free (w) ;
  1167.      
  1168.        for (i = 0; i < 5 ; i++)
  1169.          {
  1170.            printf("z%d = %+.18f %+.18f\n", i, z[2*i], z[2*i+1]) ;
  1171.          }
  1172.      }
  1173.  
  1174. The results of running the program are,
  1175.  
  1176.      bash$ ./a.out
  1177.      z0 = -0.809016994374947451 +0.587785252292473137
  1178.      z1 = -0.809016994374947451 -0.587785252292473137
  1179.      z2 = +0.309016994374947451 +0.951056516295153642
  1180.      z3 = +0.309016994374947451 -0.951056516295153642
  1181.      z4 = +1.000000000000000000 +0.000000000000000000
  1182.  
  1183. 
  1184. File: gsl-ref.info,  Node: Roots of Polynomials References and Further Reading,  Prev: Roots of Polynomials Examples,  Up: Roots of Polynomials
  1185.  
  1186. References and Further Reading
  1187. ==============================
  1188.  
  1189. The balanced-QR method and its error analysis is described in the
  1190. following papers.
  1191.  
  1192.      R.S. Martin, G. Peters and J.H. Wilkinson, "The QR Algorithm for
  1193.      Real Hessenberg Matrices", `Numerische Mathematik', 14 (1970),
  1194.      219-231.
  1195.  
  1196.      B.N. Parlett and C. Reinsch, "Balancing a Matrix for Calculation of
  1197.      Eigenvalues and Eigenvectors", `Numerische Mathematik', 13 (1969),
  1198.      293-304.
  1199.  
  1200.      A. Edelman and H. Murakami, "Polynomial roots from companion matrix
  1201.      eigenvalues", `Mathematics of Computation', Vol. 64 No. 210
  1202.      (1995), 763-776.
  1203.  
  1204.